Contracting cube

Contracting cubeΒΆ

In this demo we simulate a unit cube that is fixed at \(x = 0\) and free at \(x = 1\). We use a transversally isotropic material with fiber oriented in the \(x\)-direction.

import dolfin
try:
    from dolfin_adjoint import (
        Constant,
        DirichletBC,
        Expression,
        UnitCubeMesh,
        interpolate,
        Mesh,
    )
except ImportError:
    from dolfin import (
        Constant,
        DirichletBC,
        interpolate,
        Expression,
        UnitCubeMesh,
        Mesh,
    )
import pulse
from fenics_plotly import plot
pulse.iterate.logger.setLevel(10)
# Create mesh
N = 6
mesh = UnitCubeMesh(N, N, N)
# Create subdomains
class Free(dolfin.SubDomain):
    def inside(self, x, on_boundary):
        return x[0] > (1.0 - dolfin.DOLFIN_EPS) and on_boundary


class Fixed(dolfin.SubDomain):
    def inside(self, x, on_boundary):
        return x[0] < dolfin.DOLFIN_EPS and on_boundary
# Create a facet fuction in order to mark the subdomains
ffun = dolfin.MeshFunction("size_t", mesh, 2)
ffun.set_all(0)
# Mark the first subdomain with value 1
fixed = Fixed()
fixed_marker = 1
fixed.mark(ffun, fixed_marker)
# Mark the second subdomain with value 2
free = Free()
free_marker = 2
free.mark(ffun, free_marker)
# Create a cell function (but we are not using it)
cfun = dolfin.MeshFunction("size_t", mesh, 3)
cfun.set_all(0)
# Collect the functions containing the markers
marker_functions = pulse.MarkerFunctions(ffun=ffun, cfun=cfun)
# Create mictrotructure
V_f = dolfin.VectorFunctionSpace(mesh, "CG", 1)
# Fibers
f0 = interpolate(Expression(("1.0", "0.0", "0.0"), degree=1), V_f)
# Sheets
s0 = interpolate(Expression(("0.0", "1.0", "0.0"), degree=1), V_f)
# Fiber-sheet normal
n0 = interpolate(Expression(("0.0", "0.0", "1.0"), degree=1), V_f)
# Collect the mictrotructure
microstructure = pulse.Microstructure(f0=f0, s0=s0, n0=n0)
# Create the geometry
geometry = pulse.Geometry(
    mesh=mesh,
    marker_functions=marker_functions,
    microstructure=microstructure,
)
# Use the default material parameters
material_parameters = pulse.HolzapfelOgden.default_parameters()
# Select model for active contraction
active_model = pulse.ActiveModels.active_strain
# active_model = "active_stress"
# Set the activation
activation = Constant(0.0)
# Create material
material = pulse.HolzapfelOgden(
    active_model=active_model,
    parameters=material_parameters,
    activation=activation,
)
# Make Dirichlet boundary conditions
def dirichlet_bc(W):
    V = W if W.sub(0).num_sub_spaces() == 0 else W.sub(0)
    return DirichletBC(V, Constant((0.0, 0.0, 0.0)), fixed)
# Make Neumann boundary conditions
neumann_bc = pulse.NeumannBC(traction=Constant(0.0), marker=free_marker)
# Collect Boundary Conditions
bcs = pulse.BoundaryConditions(dirichlet=(dirichlet_bc,), neumann=(neumann_bc,))
# Create problem
problem = pulse.MechanicsProblem(geometry, material, bcs)
# Solve problem
pulse.iterate.iterate(problem, activation, 0.1)
2022-03-22 21:23:05,420 [808] DEBUG    pulse.iterate: Control: [0.0]
2022-03-22 21:23:05,420 [808] DEBUG    pulse.iterate: Target: [0.1]
2022-03-22 21:23:05,422 [808] DEBUG    pulse.iterate: Intial number of steps: 5 with step size 0.02
2022-03-22 21:23:05,423 [808] INFO     pulse.iterate: Iterating to target control (f_21)...
2022-03-22 21:23:05,423 [808] INFO     pulse.iterate: Current control: f_21 = 0.000
2022-03-22 21:23:05,424 [808] INFO     pulse.iterate: Target: 0.100
2022-03-22 21:23:05,425 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:23:05,426 [808] DEBUG    pulse.iterate: Maximum difference: 1.000e-01
2022-03-22 21:23:05,428 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:23:05,429 [808] DEBUG    pulse.iterate: Current control: 0.020
2022-03-22 21:24:52,049 [808] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2022-03-22 21:24:52,050 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:52,051 [808] DEBUG    pulse.iterate: Maximum difference: 8.000e-02
2022-03-22 21:24:52,052 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:52,053 [808] DEBUG    pulse.iterate: Maximum difference: 8.000e-02
2022-03-22 21:24:52,053 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:24:52,054 [808] DEBUG    pulse.iterate: Current control: 0.040
2022-03-22 21:24:54,609 [808] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2022-03-22 21:24:54,610 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:54,611 [808] DEBUG    pulse.iterate: Maximum difference: 6.000e-02
2022-03-22 21:24:54,612 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:54,613 [808] DEBUG    pulse.iterate: Maximum difference: 6.000e-02
2022-03-22 21:24:54,614 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:24:54,615 [808] DEBUG    pulse.iterate: Current control: 0.060
2022-03-22 21:24:55,892 [808] DEBUG    pulse.iterate: Solver did not converge...
2022-03-22 21:24:55,894 [808] DEBUG    pulse.iterate: 
NOT CONVERGING
2022-03-22 21:24:55,895 [808] DEBUG    pulse.iterate: Reduce control step
2022-03-22 21:24:55,895 [808] DEBUG    pulse.iterate: Assign old state
2022-03-22 21:24:55,901 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:55,902 [808] DEBUG    pulse.iterate: Maximum difference: 6.000e-02
2022-03-22 21:24:55,902 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:24:55,903 [808] DEBUG    pulse.iterate: Current control: 0.050
*** Warning: PETSc SNES solver diverged in 3 iterations with divergence reason DIVERGED_FNORM_NAN.
2022-03-22 21:24:57,477 [808] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2022-03-22 21:24:57,478 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:57,480 [808] DEBUG    pulse.iterate: Maximum difference: 5.000e-02
2022-03-22 21:24:57,480 [808] DEBUG    pulse.iterate: Adapt step size. New step size: 0.015
2022-03-22 21:24:57,482 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:57,483 [808] DEBUG    pulse.iterate: Maximum difference: 5.000e-02
2022-03-22 21:24:57,484 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:24:57,486 [808] DEBUG    pulse.iterate: Current control: 0.065
2022-03-22 21:24:58,769 [808] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2022-03-22 21:24:58,770 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:58,771 [808] DEBUG    pulse.iterate: Maximum difference: 3.500e-02
2022-03-22 21:24:58,772 [808] DEBUG    pulse.iterate: Adapt step size. New step size: 0.022
2022-03-22 21:24:58,774 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:58,775 [808] DEBUG    pulse.iterate: Maximum difference: 3.500e-02
2022-03-22 21:24:58,776 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:24:58,776 [808] DEBUG    pulse.iterate: Current control: 0.087
2022-03-22 21:24:59,739 [808] DEBUG    pulse.iterate: Solver did not converge...
2022-03-22 21:24:59,741 [808] DEBUG    pulse.iterate: 
NOT CONVERGING
2022-03-22 21:24:59,741 [808] DEBUG    pulse.iterate: Reduce control step
2022-03-22 21:24:59,742 [808] DEBUG    pulse.iterate: Assign old state
2022-03-22 21:24:59,748 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:24:59,749 [808] DEBUG    pulse.iterate: Maximum difference: 3.500e-02
2022-03-22 21:24:59,750 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:24:59,751 [808] DEBUG    pulse.iterate: Current control: 0.076
*** Warning: PETSc SNES solver diverged in 3 iterations with divergence reason DIVERGED_DTOL.
2022-03-22 21:25:01,011 [808] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2022-03-22 21:25:01,012 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:25:01,013 [808] DEBUG    pulse.iterate: Maximum difference: 2.375e-02
2022-03-22 21:25:01,013 [808] DEBUG    pulse.iterate: Adapt step size. New step size: 0.017
2022-03-22 21:25:01,014 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:25:01,015 [808] DEBUG    pulse.iterate: Maximum difference: 2.375e-02
2022-03-22 21:25:01,016 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:25:01,016 [808] DEBUG    pulse.iterate: Current control: 0.093
2022-03-22 21:25:02,301 [808] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2022-03-22 21:25:02,301 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:25:02,302 [808] DEBUG    pulse.iterate: Maximum difference: 6.875e-03
2022-03-22 21:25:02,303 [808] DEBUG    pulse.iterate: Adapt step size. New step size: 0.025
2022-03-22 21:25:02,304 [808] DEBUG    pulse.iterate: Check target reached: NO
2022-03-22 21:25:02,305 [808] DEBUG    pulse.iterate: Maximum difference: 6.875e-03
2022-03-22 21:25:02,305 [808] DEBUG    pulse.iterate: Change step size for final iteration
2022-03-22 21:25:02,306 [808] DEBUG    pulse.iterate: Try new control
2022-03-22 21:25:02,306 [808] DEBUG    pulse.iterate: Current control: 0.100
2022-03-22 21:25:03,268 [808] DEBUG    pulse.iterate: 
SUCCESFULL STEP:
2022-03-22 21:25:03,269 [808] DEBUG    pulse.iterate: Check target reached: YES!
2022-03-22 21:25:03,271 [808] DEBUG    pulse.iterate: Check target reached: YES!
([Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 47),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 150),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 223),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 313),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 354),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 428),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 469),
  Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 502)],
 [Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 45),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 148),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 221),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 311),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 352),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 426),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 467),
  Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 500)])
# Get displacement and hydrostatic pressure
u, p = problem.state.split(deepcopy=True)
V = dolfin.VectorFunctionSpace(mesh, "CG", 1)
u_int = interpolate(u, V)
new_mesh = Mesh(mesh)
dolfin.ALE.move(new_mesh, u_int)
fig = plot(mesh, show=False)
fig.add_plot(plot(new_mesh, color="red", show=False))
fig.show()